satkit 0.16.2

Satellite Toolkit
Documentation
{
 "cells": [
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "# Two-Line Element Set\n",
    "\n",
    "Two-Line Element sets (TLEs) are the most widely used format for distributing satellite orbital data. They encode mean orbital elements designed specifically for use with the SGP4 propagator, and are published by organizations like [CelesTrak](https://celestrak.org) and [Space-Track](https://www.space-track.org).\n",
    "\n",
    "This tutorial demonstrates:\n",
    "\n",
    "1. Loading a TLE and propagating it with SGP4 to get a position and velocity in the **TEME** (True Equator Mean Equinox) frame\n",
    "2. Rotating from TEME to the **ITRF** (Earth-fixed) frame to obtain geodetic coordinates\n",
    "3. Plotting the satellite ground track over multiple orbits"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": "## Generate State Vector\n\nWe fetch a current TLE for the International Space Station directly from CelesTrak using `sk.TLE.from_url` — no external HTTP library required. SGP4 outputs position and velocity in the TEME frame, a pseudo-inertial frame that does not account for precession or nutation. To get a location on the Earth's surface, we rotate into the ITRF frame using a quaternion from `satkit.frametransform`, then convert the Cartesian position to geodetic coordinates."
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": "import satkit as sk\n\n# Fetch the current TLE for the ISS directly from CelesTrak.  sk.TLE.from_url\n# returns a single TLE when the response contains exactly one element, or a\n# list of TLEs otherwise.  No external HTTP library is required.\nurl = \"https://celestrak.org/NORAD/elements/gp.php?CATNR=25544&FORMAT=TLE\"\niss = sk.TLE.from_url(url)\n\n# Evaluate the orbital state 6 hours after the TLE epoch\nthetime = iss.epoch + sk.duration(hours=6)\n\n# The state is output in the \"TEME\" frame, which is an approximate inertial\n# frame that does not include precession or nutation\n# pTEME is geocentric position in meters\n# vTEME is geocentric velocity in meters / second\n# for now we will ignore the velocity\npTEME, _vTEME = sk.sgp4(iss, thetime)\n\n# Suppose we want current latitude, longitude, and altitude of satellite:\n# we need to rotate into an Earth-fixed frame, the ITRF\n# We use a \"quaternion\" to represent the rotation.  Quaternion rotations\n# in the satkit toolbox can be represented as multiplications of a 3-vector\npITRF = sk.frametransform.qteme2itrf(thetime) * pTEME\n\n# Now lets make a \"ITRFCoord\" object to extract geodetic coordinates\ncoord = sk.itrfcoord(pITRF)\n\n# Get the latitude, longitude, and\n# altitude (height above ellipsoid, or hae) of the satellite\nprint(coord)"
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "## Plot Satellite Ground Track\n",
    "\n",
    "The ground track is the projection of the satellite's orbit onto the Earth's surface. The sinusoidal pattern results from the satellite's inclined orbit combined with the Earth's rotation. The `mean_motion` field in the TLE gives the number of orbits per day, which we use to compute the total time span for 5 complete orbits."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": "import satkit as sk\nimport numpy as np\nimport matplotlib.pyplot as plt\nimport scienceplots  # noqa: F401\nplt.style.use([\"science\", \"no-latex\", \"../satkit.mplstyle\"])\n%config InlineBackend.figure_formats = ['svg']\nimport warnings\nwarnings.filterwarnings(\"ignore\", \"Downloading\")\nimport cartopy.crs as ccrs\nimport cartopy.feature as cfeature\n\n# Fetch the current TLE for the ISS directly from CelesTrak\nurl = \"https://celestrak.org/NORAD/elements/gp.php?CATNR=25544&FORMAT=TLE\"\niss = sk.TLE.from_url(url)\n\n# Start the ground-track plot at the TLE epoch\nthetime = iss.epoch\n\n# plot for 5 orbits.  The mean motion in the TLE is number of orbits in a day\ntimearr = np.array([thetime + sk.duration(days=x) for x in np.linspace(0, 5/iss.mean_motion, 1000)])\n\n# Get position in the TEME frame\npTEME, _vTEME = sk.sgp4(iss, timearr)\nqarr = sk.frametransform.qteme2itrf(timearr)\npITRF = np.array([q * p for q, p in zip(qarr, pTEME)])\ncoord = [sk.itrfcoord(p) for p in pITRF]\nlat, lon, alt = zip(*[(c.latitude_deg, c.longitude_deg, c.altitude) for c in coord])\n\nfig, ax = plt.subplots(figsize=(10, 5), subplot_kw={\"projection\": ccrs.PlateCarree()})\nax.add_feature(cfeature.LAND, facecolor=\"lightgray\")\nax.add_feature(cfeature.BORDERS, linewidth=0.5)\nax.add_feature(cfeature.COASTLINE, linewidth=0.5)\n# Break line at date line crossings\nlon_arr, lat_arr = np.array(lon), np.array(lat)\nbreaks = np.where(np.abs(np.diff(lon_arr)) > 180)[0] + 1\nlon_segs = np.split(lon_arr, breaks)\nlat_segs = np.split(lat_arr, breaks)\nfor lo, la in zip(lon_segs, lat_segs):\n    ax.plot(lo, la, linewidth=1, color=\"C0\", transform=ccrs.PlateCarree())\nax.set_title(\"Ground Track\")\nax.set_global()\nplt.tight_layout()\nplt.show()\n\nfig, ax = plt.subplots(figsize=(10, 4))\nax.plot([t.as_datetime() for t in timearr], np.array(alt)/1e3)\nax.set_xlabel(\"Time\")\nax.set_ylabel(\"Altitude (km)\")\nax.set_title(\"Altitude vs Time\")\nfig.autofmt_xdate()\nplt.tight_layout()\nplt.show()"
  }
 ],
 "metadata": {
  "kernelspec": {
   "display_name": "Python 3",
   "language": "python",
   "name": "python3"
  },
  "language_info": {
   "codemirror_mode": {
    "name": "ipython",
    "version": 3
   },
   "file_extension": ".py",
   "mimetype": "text/x-python",
   "name": "python",
   "nbconvert_exporter": "python",
   "pygments_lexer": "ipython3",
   "version": "3.13.0"
  }
 },
 "nbformat": 4,
 "nbformat_minor": 2
}